Numerical Linear Algebra

Biostat/Biomath M257

Author

Dr. Hua Zhou @ UCLA

Published

April 9, 2024

System information (for reproducibility):

versioninfo()
Julia Version 1.10.2
Commit bd47eca2c8a (2024-03-01 10:14 UTC)
Build Info:
  Official https://julialang.org/ release
Platform Info:
  OS: macOS (arm64-apple-darwin22.4.0)
  CPU: 12 × Apple M2 Max
  WORD_SIZE: 64
  LIBM: libopenlibm
  LLVM: libLLVM-15.0.7 (ORCJIT, apple-m1)
Threads: 8 default, 0 interactive, 4 GC (on 8 virtual cores)
Environment:
  JULIA_NUM_THREADS = 8
  JULIA_EDITOR = code

Load packages:

using Pkg

Pkg.activate(pwd())
Pkg.instantiate()
Pkg.status()
  Activating project at `~/Documents/github.com/ucla-biostat-257/2024spring/slides/08-numalgintro`
Status `~/Documents/github.com/ucla-biostat-257/2024spring/slides/08-numalgintro/Project.toml`
  [6e4b80f9] BenchmarkTools v1.5.0
  [0e44f5e4] Hwloc v3.0.1
  [bdcacae8] LoopVectorization v0.12.169
  [6f49c342] RCall v0.14.1
  [37e2e46d] LinearAlgebra
  [9a3f8284] Random

1 Introduction

  • Topics in numerical algebra:

    • BLAS
    • solve linear equations \(\mathbf{A} \mathbf{x} = \mathbf{b}\)
    • regression computations \(\mathbf{X}^T \mathbf{X} \beta = \mathbf{X}^T \mathbf{y}\)
    • eigen-problems \(\mathbf{A} \mathbf{x} = \lambda \mathbf{x}\)
    • generalized eigen-problems \(\mathbf{A} \mathbf{x} = \lambda \mathbf{B} \mathbf{x}\)
    • singular value decompositions \(\mathbf{A} = \mathbf{U} \Sigma \mathbf{V}^T\)
    • iterative methods for numerical linear algebra
  • Except for the iterative methods, most of these numerical linear algebra tasks are implemented in the BLAS and LAPACK libraries. They form the building blocks of most statistical computing tasks (optimization, MCMC).

  • Our major goal (or learning objectives) is to

    1. know the complexity (flop count) of each task
    2. be familiar with the BLAS and LAPACK functions (what they do)
    3. do not re-invent wheels by implementing these dense linear algebra subroutines by yourself
    4. understand the need for iterative methods
    5. apply appropriate numerical algebra tools to various statistical problems
  • All high-level languages (Julia, Matlab, Python, R) call BLAS and LAPACK for numerical linear algebra.

    • Julia offers more flexibility by exposing interfaces to many BLAS/LAPACK subroutines directly. See documentation.

2 BLAS

  • BLAS stands for basic linear algebra subprograms.

  • See netlib for a complete list of standardized BLAS functions.

  • There are many implementations of BLAS.

    • Netlib provides a reference implementation.
    • Matlab uses Intel’s MKL (mathematical kernel libaries). MKL implementation is the gold standard on market. It is not open source but the compiled library is free for Linux and MacOS. However, not surprisingly, it only works on Intel CPUs.
    • Julia uses OpenBLAS. OpenBLAS is the best cross-platform, open source implementation. With the MKL.jl package, it’s also very easy to use MKL in Julia.
  • There are 3 levels of BLAS functions.

Level Example Operation Name Dimension Flops
1 \(\alpha \gets \mathbf{x}^T \mathbf{y}\) dot product \(\mathbf{x}, \mathbf{y} \in \mathbb{R}^n\) \(2n\)
1 \(\mathbf{y} \gets \mathbf{y} + \alpha \mathbf{x}\) axpy \(\alpha \in \mathbb{R}\), \(\mathbf{x}, \mathbf{y} \in \mathbb{R}^n\) \(2n\)
2 \(\mathbf{y} \gets \mathbf{y} + \mathbf{A} \mathbf{x}\) gaxpy \(\mathbf{A} \in \mathbb{R}^{m \times n}\), \(\mathbf{x} \in \mathbb{R}^n\), \(\mathbf{y} \in \mathbb{R}^m\) \(2mn\)
2 \(\mathbf{A} \gets \mathbf{A} + \mathbf{y} \mathbf{x}^T\) rank one update \(\mathbf{A} \in \mathbb{R}^{m \times n}\), \(\mathbf{x} \in \mathbb{R}^n\), \(\mathbf{y} \in \mathbb{R}^m\) \(2mn\)
3 \(\mathbf{C} \gets \mathbf{C} + \mathbf{A} \mathbf{B}\) matrix multiplication \(\mathbf{A} \in \mathbb{R}^{m \times p}\), \(\mathbf{B} \in \mathbb{R}^{p \times n}\), \(\mathbf{C} \in \mathbb{R}^{m \times n}\) \(2mnp\)
  • Typical BLAS functions support single precision (S), double precision (D), complex (C), and double complex (Z).

3 Examples

The form of a mathematical expression and the way the expression should be evaluated in actual practice may be quite different.

Some operations appear as level-3 but indeed are level-2.

Example 1. A common operation in statistics is column scaling or row scaling \[ \begin{eqnarray*} \mathbf{A} &=& \mathbf{A} \mathbf{D} \quad \text{(column scaling)} \\ \mathbf{A} &=& \mathbf{D} \mathbf{A} \quad \text{(row scaling)}, \end{eqnarray*} \] where \(\mathbf{D}\) is diagonal. For example, in generalized linear models (GLMs), the Fisher information matrix takes the form
\[ \mathbf{X}^T \mathbf{W} \mathbf{X}, \] where \(\mathbf{W}\) is a diagonal matrix with observation weights on diagonal.

Column and row scalings are essentially level-2 operations!

using BenchmarkTools, LinearAlgebra, Random

Random.seed!(257) # seed
n = 2000
A = rand(n, n) # n-by-n matrix
d = rand(n)  # n vector
D = Diagonal(d) # diagonal matrix with d as diagonal
2000×2000 Diagonal{Float64, Vector{Float64}}:
 0.0416032   ⋅         ⋅         ⋅       …   ⋅         ⋅         ⋅ 
  ⋅         0.887679   ⋅         ⋅           ⋅         ⋅         ⋅ 
  ⋅          ⋅        0.102233   ⋅           ⋅         ⋅         ⋅ 
  ⋅          ⋅         ⋅        0.08407      ⋅         ⋅         ⋅ 
  ⋅          ⋅         ⋅         ⋅           ⋅         ⋅         ⋅ 
  ⋅          ⋅         ⋅         ⋅       …   ⋅         ⋅         ⋅ 
  ⋅          ⋅         ⋅         ⋅           ⋅         ⋅         ⋅ 
  ⋅          ⋅         ⋅         ⋅           ⋅         ⋅         ⋅ 
  ⋅          ⋅         ⋅         ⋅           ⋅         ⋅         ⋅ 
  ⋅          ⋅         ⋅         ⋅           ⋅         ⋅         ⋅ 
 ⋮                                       ⋱                      
  ⋅          ⋅         ⋅         ⋅           ⋅         ⋅         ⋅ 
  ⋅          ⋅         ⋅         ⋅           ⋅         ⋅         ⋅ 
  ⋅          ⋅         ⋅         ⋅           ⋅         ⋅         ⋅ 
  ⋅          ⋅         ⋅         ⋅           ⋅         ⋅         ⋅ 
  ⋅          ⋅         ⋅         ⋅       …   ⋅         ⋅         ⋅ 
  ⋅          ⋅         ⋅         ⋅           ⋅         ⋅         ⋅ 
  ⋅          ⋅         ⋅         ⋅          0.213471   ⋅         ⋅ 
  ⋅          ⋅         ⋅         ⋅           ⋅        0.870533   ⋅ 
  ⋅          ⋅         ⋅         ⋅           ⋅         ⋅        0.318876
Dfull = diagm(d) # convert to full matrix
2000×2000 Matrix{Float64}:
 0.0416032  0.0       0.0       0.0      …  0.0       0.0       0.0
 0.0        0.887679  0.0       0.0         0.0       0.0       0.0
 0.0        0.0       0.102233  0.0         0.0       0.0       0.0
 0.0        0.0       0.0       0.08407     0.0       0.0       0.0
 0.0        0.0       0.0       0.0         0.0       0.0       0.0
 0.0        0.0       0.0       0.0      …  0.0       0.0       0.0
 0.0        0.0       0.0       0.0         0.0       0.0       0.0
 0.0        0.0       0.0       0.0         0.0       0.0       0.0
 0.0        0.0       0.0       0.0         0.0       0.0       0.0
 0.0        0.0       0.0       0.0         0.0       0.0       0.0
 ⋮                                       ⋱                      
 0.0        0.0       0.0       0.0         0.0       0.0       0.0
 0.0        0.0       0.0       0.0         0.0       0.0       0.0
 0.0        0.0       0.0       0.0         0.0       0.0       0.0
 0.0        0.0       0.0       0.0         0.0       0.0       0.0
 0.0        0.0       0.0       0.0      …  0.0       0.0       0.0
 0.0        0.0       0.0       0.0         0.0       0.0       0.0
 0.0        0.0       0.0       0.0         0.213471  0.0       0.0
 0.0        0.0       0.0       0.0         0.0       0.870533  0.0
 0.0        0.0       0.0       0.0         0.0       0.0       0.318876
# this is calling BLAS routine for matrix multiplication: O(n^3) flops
# this is SLOW!
@benchmark $A * $Dfull
BenchmarkTools.Trial: 30 samples with 1 evaluation.
 Range (minmax):  110.016 ms330.218 ms   GC (min … max): 0.00% … 0.00%
 Time  (median):     160.076 ms                GC (median):    0.00%
 Time  (mean ± σ):   169.199 ms ±  61.926 ms   GC (mean ± σ):  0.28% ± 0.57%
  █            ▁                                                 
  █▄▄▁▁▁▁▄▄▄▄▁▇█▁▁▄▁▁▄▁▄▁▁▄▁▁▁▁▁▄▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▄▁▇ ▁
  110 ms           Histogram: frequency by time          330 ms <
 Memory estimate: 30.52 MiB, allocs estimate: 2.
# dispatch to special method for diagonal matrix multiplication.
# columnwise scaling: O(n^2) flops
@benchmark $A * $D
BenchmarkTools.Trial: 1786 samples with 1 evaluation.
 Range (minmax):  2.240 ms  4.773 ms   GC (min … max):  0.00% … 32.09%
 Time  (median):     2.523 ms                GC (median):     0.00%
 Time  (mean ± σ):   2.797 ms ± 532.686 μs   GC (mean ± σ):  10.35% ± 13.58%
        ▅█▃▃                                                   
  ▂▃▃▅▅▇█████▆▄▄▃▃▂▂▂▂▂▁▁▂▁▁▂▁▂▁▁▁▁▁▁▁▁▂▂▂▂▃▃▄▅▆▆▅▄▃▃▃▃▃▂▂▂ ▃
  2.24 ms         Histogram: frequency by time        3.98 ms <
 Memory estimate: 30.52 MiB, allocs estimate: 2.
# Or we can use broadcasting (with recycling)
@benchmark $A .* transpose($d)
BenchmarkTools.Trial: 1827 samples with 1 evaluation.
 Range (minmax):  2.256 ms  4.258 ms   GC (min … max):  0.00% … 32.83%
 Time  (median):     2.460 ms                GC (median):     0.00%
 Time  (mean ± σ):   2.734 ms ± 527.168 μs   GC (mean ± σ):  10.63% ± 13.87%
    ▁█▆▆▆▂                                                    
  ▂▃███████▇▆▄▃▃▂▃▂▂▂▂▂▂▁▂▁▁▂▁▁▁▁▁▁▁▁▁▁▁▂▁▁▁▂▃▄▅▅▆▆▆▅▄▃▃▃▃▃ ▃
  2.26 ms         Histogram: frequency by time        3.83 ms <
 Memory estimate: 30.52 MiB, allocs estimate: 2.
# in-place: avoid allocate space for result
# rmul!: compute matrix-matrix product A*B, overwriting A, and return the result.
@benchmark rmul!($A, $D)
BenchmarkTools.Trial: 3314 samples with 1 evaluation.
 Range (minmax):  1.287 ms  4.918 ms   GC (min … max): 0.00% … 0.00%
 Time  (median):     1.504 ms                GC (median):    0.00%
 Time  (mean ± σ):   1.506 ms ± 130.657 μs   GC (mean ± σ):  0.00% ± 0.00%
                      ▅                                      
  ▂▁▂▂▁▂▂▂▂▃▄▄▅▅▆▅▆▆▇██▆▄▄▃▃▃▂▂▂▂▂▂▂▂▂▂▂▂▁▂▂▁▂▂▂▂▂▂▂▂▁▁▂▂▂▂ ▃
  1.29 ms         Histogram: frequency by time        1.91 ms <
 Memory estimate: 0 bytes, allocs estimate: 0.
# In-place broadcasting 
@benchmark $A .= $A .* transpose($d)
BenchmarkTools.Trial: 3051 samples with 1 evaluation.
 Range (minmax):  1.525 ms  8.799 ms   GC (min … max): 0.00% … 0.00%
 Time  (median):     1.574 ms                GC (median):    0.00%
 Time  (mean ± σ):   1.636 ms ± 326.920 μs   GC (mean ± σ):  0.00% ± 0.00%
  ▆█▃▂▁                                                    ▁
  ██████▇▆▇▆▆▄▆▅▅▆▃▅▆▃▅▅▆▄▅▃▁▄▄▄▅▅▄▄▅▄▃▃▅▅▃▃▄▃▄▁▃▃▁▅▁▃▁▁▅▄ █
  1.53 ms      Histogram: log(frequency) by time      2.99 ms <
 Memory estimate: 0 bytes, allocs estimate: 0.

Exercise: Try @turbo (SIMD) and @tturbo (SIMD) from LoopVectorization.jl package.

Note: In R or Matlab, diag(d) will create a full matrix. Be cautious using diag function: do we really need a full diagonal matrix?

using RCall

R"""
d <- runif(5)
diag(d)
"""
RObject{RealSxp}
          [,1]      [,2]      [,3]       [,4]      [,5]
[1,] 0.9431375 0.0000000 0.0000000 0.00000000 0.0000000
[2,] 0.0000000 0.1546795 0.0000000 0.00000000 0.0000000
[3,] 0.0000000 0.0000000 0.3322724 0.00000000 0.0000000
[4,] 0.0000000 0.0000000 0.0000000 0.08111416 0.0000000
[5,] 0.0000000 0.0000000 0.0000000 0.00000000 0.9322544
# This works only when Matlab is installed
using MATLAB

mat"""
d = rand(5, 1)
diag(d)
"""
ArgumentError: ArgumentError: Package MATLAB not found in current path.
- Run `import Pkg; Pkg.add("MATLAB")` to install the MATLAB package.

Example 2. Innter product between two matrices \(\mathbf{A}, \mathbf{B} \in \mathbb{R}^{m \times n}\) is often written as \[ \text{trace}(\mathbf{A}^T \mathbf{B}), \text{trace}(\mathbf{B} \mathbf{A}^T), \text{trace}(\mathbf{A} \mathbf{B}^T), \text{ or } \text{trace}(\mathbf{B}^T \mathbf{A}). \] They appear as level-3 operation (matrix multiplication with \(O(m^2n)\) or \(O(mn^2)\) flops).

Random.seed!(123)
n = 2000
A, B = randn(n, n), randn(n, n)

# slow way to evaluate tr(A'B): 2mn^2 flops
@benchmark tr(transpose($A) * $B)
BenchmarkTools.Trial: 62 samples with 1 evaluation.
 Range (minmax):  50.385 ms147.322 ms   GC (min … max): 0.00% … 0.00%
 Time  (median):     66.011 ms                GC (median):    0.00%
 Time  (mean ± σ):   81.664 ms ±  31.991 ms   GC (mean ± σ):  0.58% ± 1.40%
  █ ▁                                                           
  █▅█▁▅▇▁▁▁▅▁▅▁▅▁▁▅▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▅▁▁▁▅▅▅▁▅▅▁▁▇▅▇▇▇▇▇▇▁▇▁▇▁▅ ▁
  50.4 ms       Histogram: log(frequency) by time       130 ms <
 Memory estimate: 30.52 MiB, allocs estimate: 2.

But \(\text{trace}(\mathbf{A}^T \mathbf{B}) = <\text{vec}(\mathbf{A}), \text{vec}(\mathbf{B})>\). The latter is level-1 BLAS operation with \(O(mn)\) flops.

# smarter way to evaluate tr(A'B): 2mn flops
@benchmark dot($A, $B)
BenchmarkTools.Trial: 2698 samples with 1 evaluation.
 Range (minmax):  1.793 ms 2.206 ms   GC (min … max): 0.00% … 0.00%
 Time  (median):     1.836 ms               GC (median):    0.00%
 Time  (mean ± σ):   1.850 ms ± 48.023 μs   GC (mean ± σ):  0.00% ± 0.00%
     ▂▅███▇▇▆▆▅▅▅▃▃▂▂                                      ▂
  ▃▃▆███████████████████▇▆▆▆▇▅▅▆▅▆▅▅▅▇▅▅▆▆▅▆▅▆▆█▆▆▆▅▆▇▅▆▄▆ █
  1.79 ms      Histogram: log(frequency) by time     2.05 ms <
 Memory estimate: 0 bytes, allocs estimate: 0.

Example 3. Similarly \(\text{diag}(\mathbf{A}^T \mathbf{B})\) can be calculated in \(O(mn)\) flops.

# slow way to evaluate diag(A'B): O(n^3)
@benchmark diag(transpose($A) * $B)
BenchmarkTools.Trial: 93 samples with 1 evaluation.
 Range (minmax):  50.588 ms75.902 ms   GC (min … max): 0.00% … 0.00%
 Time  (median):     52.324 ms               GC (median):    0.00%
 Time  (mean ± σ):   53.967 ms ±  4.768 ms   GC (mean ± σ):  0.89% ± 1.54%
  █ ▁                                                      
  █▇█▇█▃▃▄▃▁▁▃▃▃▃▄▃▃▃▁▃▃▁▁▁▁▁▃▁▁▃▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▃ ▁
  50.6 ms         Histogram: frequency by time          74 ms <
 Memory estimate: 30.53 MiB, allocs estimate: 3.
# smarter way to evaluate diag(A'B): O(n^2)
@benchmark Diagonal(vec(sum($A .* $B, dims = 1)))
BenchmarkTools.Trial: 1257 samples with 1 evaluation.
 Range (minmax):  3.330 ms 15.680 ms   GC (min … max): 0.00% …  9.60%
 Time  (median):     3.640 ms                GC (median):    0.00%
 Time  (mean ± σ):   3.972 ms ± 774.977 μs   GC (mean ± σ):  9.38% ± 12.67%
   ▄█▅▃▁                                                       
  ▅██████▇▇▅▅▄▃▃▂▂▂▂▂▂▁▂▂▁▁▂▁▂▃▅▆▃▄▄▄▄▄▃▃▂▃▃▂▂▂▂▂▂▂▂▂▃▃▃▃▂▂ ▃
  3.33 ms         Histogram: frequency by time        5.93 ms <
 Memory estimate: 30.53 MiB, allocs estimate: 5.

To get rid of allocation of intermediate arrays at all, we can just write a double loop or use dot function.

function diag_matmul!(d, A, B)
    m, n = size(A)
    @assert size(B) == (m, n) "A and B should have same size"
    fill!(d, 0)
    for j in 1:n, i in 1:m
        d[j] += A[i, j] * B[i, j]
    end
    Diagonal(d)
end

d = zeros(eltype(A), size(A, 2))
@benchmark diag_matmul!($d, $A, $B)
BenchmarkTools.Trial: 1465 samples with 1 evaluation.
 Range (minmax):  3.319 ms 13.782 ms   GC (min … max): 0.00% … 0.00%
 Time  (median):     3.374 ms                GC (median):    0.00%
 Time  (mean ± σ):   3.409 ms ± 509.119 μs   GC (mean ± σ):  0.00% ± 0.00%
  ▅▄ ▁▃▂▃▇▇▇▇█▆▃▆▄▄ ▁                                        
  ██▇██████████████▇█▆▆▆▄▅▅▅▅▅▄▃▃▃▃▃▃▂▃▃▂▁▂▂▁▁▂▁▁▁▁▂▂▁▁▂▁▂▂ ▄
  3.32 ms         Histogram: frequency by time        3.58 ms <
 Memory estimate: 0 bytes, allocs estimate: 0.

Exercise: Try @turbo (SIMD) and @tturbo (SIMD) from LoopVectorization.jl package.

4 Memory hierarchy and level-3 fraction

Key to high performance is effective use of memory hierarchy. True on all architectures.

  • Flop count is not the sole determinant of algorithm efficiency. Another important factor is data movement through the memory hierarchy.

Source: https://cs.brown.edu/courses/csci1310/2020/assign/labs/lab4.html

  • In Julia, we can query the CPU topology by the Hwloc.jl package. For example, this laptop runs an Apple M2 Max chip with 4 efficiency cores and 8 performance cores.
using Hwloc

topology_graphical()
/------------------------------------------------------------------------------------------------------------------------------------------------------------------\
| Machine (2208MB total)                                                                                                                                           |
|                                                                                                                                                                  |
| /----------------------------------------------------------------------------------------------------------------------------------------\  /------------------\ |
| | Package L#0                                                                                                                            |  | CoProc opencl0d0 | |
| |                                                                                                                                        |  |                  | |
| | /------------------------------------------------------------------------------------------------------------------------------------\ |  | 38 compute units | |
| | | NUMANode L#0 P#0 (2208MB)                                                                                                          | |  |                  | |
| | \------------------------------------------------------------------------------------------------------------------------------------/ |  | 72 GB            | |
| |                                                                                                                                        |  \------------------/ |
| | /----------------------------------------------------------------\  /----------------------------------------------------------------\ |                       |
| | | L2 (4096KB)                                                    |  | L2 (16MB)                                                      | |                       |
| | \----------------------------------------------------------------/  \----------------------------------------------------------------/ |                       |
| |                                                                                                                                        |                       |
| | /-------------\  /-------------\  /-------------\  /-------------\  /-------------\  /-------------\  /-------------\  /-------------\ |                       |
| | | L1d (64KB)  |  | L1d (64KB)  |  | L1d (64KB)  |  | L1d (64KB)  |  | L1d (128KB) |  | L1d (128KB) |  | L1d (128KB) |  | L1d (128KB) | |                       |
| | \-------------/  \-------------/  \-------------/  \-------------/  \-------------/  \-------------/  \-------------/  \-------------/ |                       |
| |                                                                                                                                        |                       |
| | /-------------\  /-------------\  /-------------\  /-------------\  /-------------\  /-------------\  /-------------\  /-------------\ |                       |
| | | L1i (128KB) |  | L1i (128KB) |  | L1i (128KB) |  | L1i (128KB) |  | L1i (192KB) |  | L1i (192KB) |  | L1i (192KB) |  | L1i (192KB) | |                       |
| | \-------------/  \-------------/  \-------------/  \-------------/  \-------------/  \-------------/  \-------------/  \-------------/ |                       |
| |                                                                                                                                        |                       |
| | /-------------\  /-------------\  /-------------\  /-------------\  /-------------\  /-------------\  /-------------\  /-------------\ |                       |
| | | Core L#0    |  | Core L#1    |  | Core L#2    |  | Core L#3    |  | Core L#4    |  | Core L#5    |  | Core L#6    |  | Core L#7    | |                       |
| | |             |  |             |  |             |  |             |  |             |  |             |  |             |  |             | |                       |
| | | /---------\ |  | /---------\ |  | /---------\ |  | /---------\ |  | /---------\ |  | /---------\ |  | /---------\ |  | /---------\ | |                       |
| | | | PU L#0  | |  | | PU L#1  | |  | | PU L#2  | |  | | PU L#3  | |  | | PU L#4  | |  | | PU L#5  | |  | | PU L#6  | |  | | PU L#7  | | |                       |
| | | |         | |  | |         | |  | |         | |  | |         | |  | |         | |  | |         | |  | |         | |  | |         | | |                       |
| | | |   P#0   | |  | |   P#1   | |  | |   P#2   | |  | |   P#3   | |  | |   P#4   | |  | |   P#5   | |  | |   P#6   | |  | |   P#7   | | |                       |
| | | \---------/ |  | \---------/ |  | \---------/ |  | \---------/ |  | \---------/ |  | \---------/ |  | \---------/ |  | \---------/ | |                       |
| | \-------------/  \-------------/  \-------------/  \-------------/  \-------------/  \-------------/  \-------------/  \-------------/ |                       |
| |                                                                                                                                        |                       |
| | /----------------------------------------------------------------\                                                                     |                       |
| | | L2 (16MB)                                                      |                                                                     |                       |
| | \----------------------------------------------------------------/                                                                     |                       |
| |                                                                                                                                        |                       |
| | /-------------\  /-------------\  /-------------\  /-------------\                                                                     |                       |
| | | L1d (128KB) |  | L1d (128KB) |  | L1d (128KB) |  | L1d (128KB) |                                                                     |                       |
| | \-------------/  \-------------/  \-------------/  \-------------/                                                                     |                       |
| |                                                                                                                                        |                       |
| | /-------------\  /-------------\  /-------------\  /-------------\                                                                     |                       |
| | | L1i (192KB) |  | L1i (192KB) |  | L1i (192KB) |  | L1i (192KB) |                                                                     |                       |
| | \-------------/  \-------------/  \-------------/  \-------------/                                                                     |                       |
| |                                                                                                                                        |                       |
| | /-------------\  /-------------\  /-------------\  /-------------\                                                                     |                       |
| | | Core L#8    |  | Core L#9    |  | Core L#10   |  | Core L#11   |                                                                     |                       |
| | |             |  |             |  |             |  |             |                                                                     |                       |
| | | /---------\ |  | /---------\ |  | /---------\ |  | /---------\ |                                                                     |                       |
| | | | PU L#8  | |  | | PU L#9  | |  | | PU L#10 | |  | | PU L#11 | |                                                                     |                       |
| | | |         | |  | |         | |  | |         | |  | |         | |                                                                     |                       |
| | | |   P#8   | |  | |   P#9   | |  | |  P#10   | |  | |  P#11   | |                                                                     |                       |
| | | \---------/ |  | \---------/ |  | \---------/ |  | \---------/ |                                                                     |                       |
| | \-------------/  \-------------/  \-------------/  \-------------/                                                                     |                       |
| \----------------------------------------------------------------------------------------------------------------------------------------/                       |
\------------------------------------------------------------------------------------------------------------------------------------------------------------------/
  • For example, Xeon X5650 CPU has a theoretical throughput of 128 DP GFLOPS but a max memory bandwidth of 32GB/s.

  • Can we keep CPU cores busy with enough deliveries of matrix data and ship the results to memory fast enough to avoid backlog?
    Answer: use high-level BLAS as much as possible.

BLAS Dimension Mem. Refs. Flops Ratio
Level 1: \(\mathbf{y} \gets \mathbf{y} + \alpha \mathbf{x}\) \(\mathbf{x}, \mathbf{y} \in \mathbb{R}^n\) \(3n\) \(2n\) 3:2
Level 2: \(\mathbf{y} \gets \mathbf{y} + \mathbf{A} \mathbf{x}\) \(\mathbf{x}, \mathbf{y} \in \mathbb{R}^n\), \(\mathbf{A} \in \mathbb{R}^{n \times n}\) \(n^2\) \(2n^2\) 1:2
Level 3: \(\mathbf{C} \gets \mathbf{C} + \mathbf{A} \mathbf{B}\) \(\mathbf{A}, \mathbf{B}, \mathbf{C} \in\mathbb{R}^{n \times n}\) \(4n^2\) \(2n^3\) 2:n
  • Higher level BLAS (3 or 2) make more effective use of arithmetic logic units (ALU) by keeping them busy. Surface-to-volume effect.

Source: Jack Dongarra’s slides.

  • A distinction between LAPACK and LINPACK (older version of R uses LINPACK) is that LAPACK makes use of higher level BLAS as much as possible (usually by smart partitioning) to increase the so-called level-3 fraction.

  • To appreciate the efforts in an optimized BLAS implementation such as OpenBLAS (evolved from GotoBLAS), see the Quora question, especially the video. Bottomline is

Get familiar with (good implementations of) BLAS/LAPACK and use them as much as possible.

5 Effect of data layout

  • Data layout in memory affects algorithmic efficiency too. It is much faster to move chunks of data in memory than retrieving/writing scattered data.

  • Storage mode: column-major (Fortran, Matlab, R, Julia) vs row-major (C/C++).

  • Cache line is the minimum amount of cache which can be loaded and stored to memory.

    • x86 CPUs: 64 bytes
    • ARM CPUs: 32 bytes

  • In Julia, we can query the cache line size by Hwloc.jl.
# Apple Silicon (M1/M2 chips) don't have L3 cache
Hwloc.cachelinesize()
ErrorException: Your system doesn't seem to have an L3 cache.
  • Accessing column-major stored matrix by rows (\(ij\) looping) causes lots of cache misses.

  • Take matrix multiplication as an example \[ \mathbf{C} \gets \mathbf{C} + \mathbf{A} \mathbf{B}, \quad \mathbf{A} \in \mathbb{R}^{m \times p}, \mathbf{B} \in \mathbb{R}^{p \times n}, \mathbf{C} \in \mathbb{R}^{m \times n}. \] Assume the storage is column-major, such as in Julia. There are 6 variants of the algorithms according to the order in the triple loops.

    • jki or kji looping:
# inner most loop
for i in 1:m
    C[i, j] = C[i, j] + A[i, k] * B[k, j]
end
- `ikj` or `kij` looping:
# inner most loop        
for j in 1:n
    C[i, j] = C[i, j] + A[i, k] * B[k, j]
end
  • ijk or jik looping:
# inner most loop        
for k in 1:p
    C[i, j] = C[i, j] + A[i, k] * B[k, j]
end
  • We pay attention to the innermost loop, where the vector calculation occurs. The associated stride when accessing the three matrices in memory (assuming column-major storage) is
Variant A Stride B Stride C Stride
\(jki\) or \(kji\) Unit 0 Unit
\(ikj\) or \(kij\) 0 Non-Unit Non-Unit
\(ijk\) or \(jik\) Non-Unit Unit 0

Apparently the variants \(jki\) or \(kji\) are preferred.

"""
    matmul_by_loop!(A, B, C, order)

Overwrite `C` by `A * B`. `order` indicates the looping order for triple loop.
"""
function matmul_by_loop!(A::Matrix, B::Matrix, C::Matrix, order::String)
    
    m = size(A, 1)
    p = size(A, 2)
    n = size(B, 2)
    fill!(C, 0)
    
    if order == "jki"
        @inbounds for j = 1:n, k = 1:p, i = 1:m
            C[i, j] += A[i, k] * B[k, j]
        end
    end

    if order == "kji"
        @inbounds for k = 1:p, j = 1:n, i = 1:m
            C[i, j] += A[i, k] * B[k, j]
        end
    end
    
    if order == "ikj"
        @inbounds for i = 1:m, k = 1:p, j = 1:n
            C[i, j] += A[i, k] * B[k, j]
        end
    end

    if order == "kij"
        @inbounds for k = 1:p, i = 1:m, j = 1:n
            C[i, j] += A[i, k] * B[k, j]
        end
    end
    
    if order == "ijk"
        @inbounds for i = 1:m, j = 1:n, k = 1:p
            C[i, j] += A[i, k] * B[k, j]
        end
    end
    
    if order == "jik"
        @inbounds for j = 1:n, i = 1:m, k = 1:p
            C[i, j] += A[i, k] * B[k, j]
        end
    end
    
end

using Random

Random.seed!(123)
m, p, n = 2000, 100, 2000
A = rand(m, p)
B = rand(p, n)
C = zeros(m, n);
  • \(jki\) and \(kji\) looping:
using BenchmarkTools

@benchmark matmul_by_loop!($A, $B, $C, "jki")
BenchmarkTools.Trial: 86 samples with 1 evaluation.
 Range (minmax):  57.729 ms70.826 ms   GC (min … max): 0.00% … 0.00%
 Time  (median):     58.222 ms               GC (median):    0.00%
 Time  (mean ± σ):   58.433 ms ±  1.412 ms   GC (mean ± σ):  0.00% ± 0.00%
   ▁  █  ▆  ▁▁▁▁   ▁   ▁  ▆    ▁                        
  ▇█▇▁█▇▄█▇▇████▄▇▁█▄▇█▄▇▄▄▇█▁▁█▄▄▁▁█▄▄▄▄▁▇▁▁▁▁▁▄▄▁▄▁▁▁▇▄▁▄ ▁
  57.7 ms         Histogram: frequency by time        59.3 ms <
 Memory estimate: 0 bytes, allocs estimate: 0.
@benchmark matmul_by_loop!($A, $B, $C, "kji")
BenchmarkTools.Trial: 27 samples with 1 evaluation.
 Range (minmax):  183.530 ms212.516 ms   GC (min … max): 0.00% … 0.00%
 Time  (median):     184.559 ms                GC (median):    0.00%
 Time  (mean ± σ):   186.442 ms ±   5.699 ms   GC (mean ± σ):  0.00% ± 0.00%
  █▆                                                            
  ██▁▁▆▆▁▁▁▄▁▁▁▄▁▁▁▄▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▄ ▁
  184 ms           Histogram: frequency by time          213 ms <
 Memory estimate: 0 bytes, allocs estimate: 0.
  • \(ikj\) and \(kij\) looping:
@benchmark matmul_by_loop!($A, $B, $C, "ikj")
BenchmarkTools.Trial: 10 samples with 1 evaluation.
 Range (minmax):  509.252 ms527.728 ms   GC (min … max): 0.00% … 0.00%
 Time  (median):     515.239 ms                GC (median):    0.00%
 Time  (mean ± σ):   515.454 ms ±   5.554 ms   GC (mean ± σ):  0.00% ± 0.00%
  █     ▁   ▁     ▁     ▁  ▁    ▁▁                            ▁  
  █▁▁▁▁▁█▁▁▁█▁▁▁▁▁█▁▁▁█▁▁█▁▁▁▁██▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁█ ▁
  509 ms           Histogram: frequency by time          528 ms <
 Memory estimate: 0 bytes, allocs estimate: 0.
@benchmark matmul_by_loop!($A, $B, $C, "kij")
BenchmarkTools.Trial: 10 samples with 1 evaluation.
 Range (minmax):  507.229 ms530.723 ms   GC (min … max): 0.00% … 0.00%
 Time  (median):     511.375 ms                GC (median):    0.00%
 Time  (mean ± σ):   513.003 ms ±   6.944 ms   GC (mean ± σ):  0.00% ± 0.00%
  ████     █   █  ██                                      █  
  ████▁▁▁▁▁██▁▁▁▁█▁▁██▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁█ ▁
  507 ms           Histogram: frequency by time          531 ms <
 Memory estimate: 0 bytes, allocs estimate: 0.
  • \(ijk\) and \(jik\) looping:
@benchmark matmul_by_loop!($A, $B, $C, "ijk")
BenchmarkTools.Trial: 21 samples with 1 evaluation.
 Range (minmax):  244.667 ms265.676 ms   GC (min … max): 0.00% … 0.00%
 Time  (median):     247.189 ms                GC (median):    0.00%
 Time  (mean ± σ):   249.187 ms ±   5.163 ms   GC (mean ± σ):  0.00% ± 0.00%
   █▃▃  ▃                                                   
  ▇███▁▁█▇█▁▁▁▁▇▁▁▁▁▁▁▇▁▇▁▇▁▇▁▁▁▁▇▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▇ ▁
  245 ms           Histogram: frequency by time          266 ms <
 Memory estimate: 0 bytes, allocs estimate: 0.
@benchmark matmul_by_loop!($A, $B, $C, "ijk")
BenchmarkTools.Trial: 21 samples with 1 evaluation.
 Range (minmax):  245.461 ms255.071 ms   GC (min … max): 0.00% … 0.00%
 Time  (median):     247.440 ms                GC (median):    0.00%
 Time  (mean ± σ):   248.889 ms ±   3.108 ms   GC (mean ± σ):  0.00% ± 0.00%
  █▁▁   ▁▁▁▁█ ▁▁         ▁  ▁    ▁     ▁       ▁ ▁  ▁▁        ▁  
  ███▁▁▁███████▁▁▁▁▁▁▁▁█▁▁█▁▁▁▁█▁▁▁▁▁█▁▁▁▁▁▁▁█▁█▁▁██▁▁▁▁▁▁▁▁█ ▁
  245 ms           Histogram: frequency by time          255 ms <
 Memory estimate: 0 bytes, allocs estimate: 0.
  • Question: Can our loop beat BLAS? Julia wraps BLAS library for matrix multiplication. We see BLAS library wins hands down (multi-threading, Strassen algorithm, higher level-3 fraction by block outer product).
@benchmark mul!($C, $A, $B)
BenchmarkTools.Trial: 1708 samples with 1 evaluation.
 Range (minmax):  2.526 ms24.014 ms   GC (min … max): 0.00% … 0.00%
 Time  (median):     2.643 ms               GC (median):    0.00%
 Time  (mean ± σ):   2.926 ms ±  1.033 ms   GC (mean ± σ):  0.00% ± 0.00%
  ██▂▃▂▃▃▁▁▁ ▁                                             ▁
  ██████████████▇▇█▇▇▇▆▆▅▆▄▄▆▆▁▁▆▄▄▅▁▁▄▄▄▁▁▄▄▁▁▁▄▁▅▁▁▁▄▄▄▄ █
  2.53 ms      Histogram: log(frequency) by time     6.75 ms <
 Memory estimate: 0 bytes, allocs estimate: 0.
# direct call of BLAS wrapper function
@benchmark LinearAlgebra.BLAS.gemm!('N', 'N', 1.0, $A, $B, 0.0, $C)
BenchmarkTools.Trial: 1653 samples with 1 evaluation.
 Range (minmax):  2.554 ms26.166 ms   GC (min … max): 0.00% … 0.00%
 Time  (median):     2.672 ms               GC (median):    0.00%
 Time  (mean ± σ):   3.022 ms ±  1.278 ms   GC (mean ± σ):  0.00% ± 0.00%
  █▄▂▂▃▂▂▂▁  ▁                                              
  ██████████████▅▇▆▆▆▆▅▇▆▇▄▄▄▅▄▄▁▄▄▄▄▄▁▄▄▄▄▄▁▄▁▄▁▄▄▄▁▄▄▄▁▄ █
  2.55 ms      Histogram: log(frequency) by time     7.47 ms <
 Memory estimate: 0 bytes, allocs estimate: 0.

Question (again): Can our loop beat BLAS?

Exercise: Annotate the loop in matmul_by_loop! by @turbo and @tturbo (multi-threading) and benchmark again.

6 BLAS in R

  • Tip for R users. Standard R distribution from CRAN uses a very out-dated BLAS/LAPACK library.
using RCall

R"""
sessionInfo()
"""
RObject{VecSxp}
R version 4.3.2 (2023-10-31)
Platform: aarch64-apple-darwin20 (64-bit)
Running under: macOS Sonoma 14.4.1

Matrix products: default
BLAS:   /System/Library/Frameworks/Accelerate.framework/Versions/A/Frameworks/vecLib.framework/Versions/A/libBLAS.dylib 
LAPACK: /Library/Frameworks/R.framework/Versions/4.3-arm64/Resources/lib/libRlapack.dylib;  LAPACK version 3.11.0

locale:
[1] C

time zone: America/Los_Angeles
tzcode source: internal

attached base packages:
[1] stats     graphics  grDevices utils     datasets  methods   base     

loaded via a namespace (and not attached):
[1] compiler_4.3.2
R"""
library(dplyr)
library(bench)
A <- $A
B <- $B
bench::mark(A %*% B) %>%
  print(width = Inf)
""";
┌ Warning: RCall.jl: 
│ Attaching package: 'dplyr'
│ 
│ The following objects are masked from 'package:stats':
│ 
│     filter, lag
│ 
│ The following objects are masked from 'package:base':
│ 
│     intersect, setdiff, setequal, union
│ 
└ @ RCall /Users/huazhou/.julia/packages/RCall/dDAVd/src/io.jl:172
┌ Warning: RCall.jl: Warning: Some expressions had a GC in every iteration; so filtering is disabled.
└ @ RCall /Users/huazhou/.julia/packages/RCall/dDAVd/src/io.jl:172
# A tibble: 1 x 13
  expression      min   median `itr/sec` mem_alloc `gc/sec` n_itr  n_gc
  <bch:expr> <bch:tm> <bch:tm>     <dbl> <bch:byt>    <dbl> <int> <dbl>
1 A %*% B       124ms    125ms      7.75    30.5MB     7.75     4     4
  total_time result                memory             time          
    <bch:tm> <list>                <list>             <list>        
1      516ms <dbl [2,000 x 2,000]> <Rprofmem [1 x 3]> <bench_tm [4]>
  gc              
  <list>          
1 <tibble [4 x 3]>
  • Re-build R from source using OpenBLAS or MKL will immediately boost linear algebra performance in R. Google build R using MKL to get started. Similarly we can build Julia using MKL.

  • Matlab uses MKL. Usually it’s very hard to beat Matlab in terms of linear algebra.

using MATLAB

mat"""
f = @() $A * $B;
timeit(f)
"""
ArgumentError: ArgumentError: Package MATLAB not found in current path.
- Run `import Pkg; Pkg.add("MATLAB")` to install the MATLAB package.

7 Avoid memory allocation: some examples

7.1 Transposing matrix is an expensive memory operation

In R, the command

t(A) %*% x

will first transpose A then perform matrix multiplication, causing unnecessary memory allocation

using Random, LinearAlgebra, BenchmarkTools
Random.seed!(123)

n = 1000
A = rand(n, n)
x = rand(n);
R"""
A <- $A
x <- $x
bench::mark(t(A) %*% x) %>%
  print(width = Inf)
""";
# A tibble: 1 x 13
  expression      min   median `itr/sec` mem_alloc `gc/sec` n_itr  n_gc
  <bch:expr> <bch:tm> <bch:tm>     <dbl> <bch:byt>    <dbl> <int> <dbl>
1 t(A) %*% x    2.1ms   2.22ms      446.    7.64MB     91.4   166    34
  total_time result            memory             time            
    <bch:tm> <list>            <list>             <list>          
1      372ms <dbl [1,000 x 1]> <Rprofmem [2 x 3]> <bench_tm [200]>
  gc                
  <list>            
1 <tibble [200 x 3]>

Julia is avoids transposing matrix whenever possible. The Transpose type only creates a view of the transpose of matrix data.

typeof(transpose(A))
Transpose{Float64, Matrix{Float64}}
fieldnames(typeof(transpose(A)))
(:parent,)
# same data in tranpose(A) and original matrix A
pointer(transpose(A).parent), pointer(A)
(Ptr{Float64} @0x0000000129840000, Ptr{Float64} @0x0000000129840000)
# dispatch to BLAS
# does *not* actually transpose the matrix
@benchmark transpose($A) * $x
BenchmarkTools.Trial: 10000 samples with 1 evaluation.
 Range (minmax):  22.625 μs 1.402 ms   GC (min … max): 0.00% … 0.00%
 Time  (median):     22.959 μs               GC (median):    0.00%
 Time  (mean ± σ):   25.961 μs ± 24.431 μs   GC (mean ± σ):  0.00% ± 0.00%
  ▇▃                                  ▁                      ▁
  ███▇▅▆▅▅▆▇▇▇▇▇▆▆▆▆▆▆▆▆▅▅▆▅▆▅▆▅▅▅▁▅███▇▅▅▆▆▆▅▅▆▅▆▅▅▅▄▄▄▃▄▅ █
  22.6 μs      Histogram: log(frequency) by time      65.8 μs <
 Memory estimate: 8.00 KiB, allocs estimate: 1.
# pre-allocate result
out = zeros(size(A, 2))
@benchmark mul!($out, transpose($A), $x)
BenchmarkTools.Trial: 10000 samples with 1 evaluation.
 Range (minmax):  22.542 μs 2.993 ms   GC (min … max): 0.00% … 0.00%
 Time  (median):     22.834 μs               GC (median):    0.00%
 Time  (mean ± σ):   26.407 μs ± 54.119 μs   GC (mean ± σ):  0.00% ± 0.00%
  █                                         ▁           ▁
  ███▆██▆▇▇▆▇▇▇▇▆▆▅▄▄▅▄▅▄▄▄▅▅▅▄▅▄▄▄▄▁▄▁▄▅▇▄▆▆▇▇███▆▄▄▃▄▄▄▅▅ █
  22.5 μs      Histogram: log(frequency) by time      56.7 μs <
 Memory estimate: 0 bytes, allocs estimate: 0.
# or call BLAS wrapper directly
@benchmark LinearAlgebra.BLAS.gemv!('T', 1.0, $A, $x, 0.0, $out)
BenchmarkTools.Trial: 10000 samples with 1 evaluation.
 Range (minmax):  22.666 μs86.333 μs   GC (min … max): 0.00% … 0.00%
 Time  (median):     22.875 μs               GC (median):    0.00%
 Time  (mean ± σ):   23.589 μs ±  2.546 μs   GC (mean ± σ):  0.00% ± 0.00%
  █▃▃▄▂▁                                                    ▁
  ██████▆▆▅▄▁██▆█▇███▆▇▇▇█▆▇▅▅▅▆▆▆▅▆▅▅▆▅▆▆▆▆▇▅▅▅▄▄▁▃▄▄▄▅▄▅▅ █
  22.7 μs      Histogram: log(frequency) by time      35.2 μs <
 Memory estimate: 0 bytes, allocs estimate: 0.

7.2 Broadcast (dot operation) fuses loops and avoids memory allocation

Broadcasting in Julia achieves vectorized code without creating intermediate arrays.

Suppose we want to calculate elementsize maximum of absolute values of two large arrays. In R or Matlab, the command

max(abs(X), abs(Y))

will create two intermediate arrays and then one result array.

using RCall

Random.seed!(123)
X, Y = rand(1000, 1000), rand(1000, 1000)

R"""
library(dplyr)
library(bench)
bench::mark(max(abs($X), abs($Y))) %>%
  print(width = Inf)
""";
# A tibble: 1 x 13
  expression                           min   median `itr/sec` mem_alloc `gc/sec`
  <bch:expr>                      <bch:tm> <bch:tm>     <dbl> <bch:byt>    <dbl>
1 max(abs(`#JL`$X), abs(`#JL`$Y))   2.85ms   3.74ms      269.    15.3MB     227.
  n_itr  n_gc total_time result    memory             time            
  <int> <dbl>   <bch:tm> <list>    <list>             <list>          
1    63    53      234ms <dbl [1]> <Rprofmem [2 x 3]> <bench_tm [116]>
  gc                
  <list>            
1 <tibble [116 x 3]>

In Julia, dot operations are fused so no intermediate arrays are created.

# no intermediate arrays created, only result array created
@benchmark max.(abs.($X), abs.($Y))
BenchmarkTools.Trial: 5994 samples with 1 evaluation.
 Range (minmax):  230.583 μs  6.933 ms   GC (min … max):  0.00% … 55.37%
 Time  (median):     736.792 μs                GC (median):     0.00%
 Time  (mean ± σ):   832.411 μs ± 522.295 μs   GC (mean ± σ):  15.75% ± 18.35%
  ▅▅▂      ▆█▇▅▃▂                             ▁   ▁           ▂
  ███▇▆▁▃▁▅███████▆▆▃▄▄▆▇▆▆▆▅▆▇▆▆▆▅▅▅▅▅▆████▇▇██▇▇███▇█▆▆▆▆▄▅ █
  231 μs        Histogram: log(frequency) by time       2.92 ms <
 Memory estimate: 7.63 MiB, allocs estimate: 2.

Pre-allocating result array gets rid of memory allocation at all.

# no memory allocation at all!
Z = zeros(size(X)) # zero matrix of same size as X
@benchmark $Z .= max.(abs.($X), abs.($Y)) # .= (vs =) is important!
BenchmarkTools.Trial: 10000 samples with 1 evaluation.
 Range (minmax):  213.292 μs573.292 μs   GC (min … max): 0.00% … 0.00%
 Time  (median):     229.792 μs                GC (median):    0.00%
 Time  (mean ± σ):   233.820 μs ±  11.851 μs   GC (mean ± σ):  0.00% ± 0.00%
                  ▃█▇▂                                           
  ▁▃▅▄▃▂▃▃▂▁▂▂▂▂▃▆████▇▆▄▃▂▂▂▂▂▂▂▂▂▂▄█▇▄▄▅▄▃▂▂▂▂▂▂▂▁▂▁▁▁▁▁▁▁▁ ▂
  213 μs           Histogram: frequency by time          264 μs <
 Memory estimate: 0 bytes, allocs estimate: 0.

Exercise: Annotate the broadcasting by @turbo and @tturbo (multi-threading) and benchmark again.

7.3 Views

View avoids creating extra copy of matrix data.

Random.seed!(123)
A = randn(1000, 1000)

# sum entries in a sub-matrix
@benchmark sum($A[1:2:500, 1:2:500])
BenchmarkTools.Trial: 10000 samples with 1 evaluation.
 Range (minmax):  45.291 μs  4.860 ms   GC (min … max):  0.00% … 97.70%
 Time  (median):     68.750 μs                GC (median):     0.00%
 Time  (mean ± σ):   74.404 μs ± 121.537 μs   GC (mean ± σ):  10.36% ±  6.73%
  █     ▁                 ▆▃                                    
  █▅▁▁▁▂█▄▁▁▂▂▁▁▁▁▁▁▁▁▁▁▁▃██▅▅▆▇▆▇▆▆▅▅▅▄▄▄▄▃▃▂▂▂▂▂▂▂▂▂▁▁▁▁▁▁ ▂
  45.3 μs         Histogram: frequency by time         93.4 μs <
 Memory estimate: 488.42 KiB, allocs estimate: 2.
# view avoids creating a separate sub-matrix
# unfortunately it's much slower possibly because of cache misses
@benchmark sum(@view $A[1:2:500, 1:2:500])
BenchmarkTools.Trial: 10000 samples with 1 evaluation.
 Range (minmax):  54.958 μs207.333 μs   GC (min … max): 0.00% … 0.00%
 Time  (median):     55.167 μs                GC (median):    0.00%
 Time  (mean ± σ):   56.099 μs ±   3.432 μs   GC (mean ± σ):  0.00% ± 0.00%
  █▅▃▃▅▁▁     ▁▁                                             ▁
  █████████▇███████████▇▇▇▇▇▇▇▇▇▆▇▇█▇▇▆▆▆▆▆▇▆▅▆▆▅▅▅▅▅▆▅▃▄▅▅▅ █
  55 μs         Histogram: log(frequency) by time      65.8 μs <
 Memory estimate: 0 bytes, allocs estimate: 0.

The @views macro, which can be useful in some operations.

@benchmark @views sum($A[1:2:500, 1:2:500])
BenchmarkTools.Trial: 10000 samples with 1 evaluation.
 Range (minmax):  54.958 μs 4.078 ms   GC (min … max): 0.00% … 0.00%
 Time  (median):     55.250 μs               GC (median):    0.00%
 Time  (mean ± σ):   58.562 μs ± 69.852 μs   GC (mean ± σ):  0.00% ± 0.00%
  █▅▃▄▂    ▁▁▁ ▂▁  ▁                                        ▁
  █████▇▆▇▆███▇██████▇▆▇▆▆▆▇▆▆▆▇▇▇▇▆▆▅▆▇▅▆▆▅▄▄▅▄▄▅▅▅▃▅▅▄▄▃▅ █
  55 μs        Histogram: log(frequency) by time        69 μs <
 Memory estimate: 0 bytes, allocs estimate: 0.

Exercise: Here we saw that, although we avoids memory allocation, the actual run times are longer. Why?